Example for cross sectional flux method¶
Author: Gerrit Kuhlmann (gerrit.kuhlmann@empa.ch)
This example code applies the cross sectional flux method to TROPOMI NO2 observations to quantify the emissions of the Janschwalde power plant.
%matplotlib inline
import os
import numpy as np
import pandas as pd
import scipy.stats
import xarray as xr
import ddeq
import ucat
Create source dataset and setup domain for plotting¶
# Source dataset
sources = xr.Dataset({
"source": xr.DataArray(["6298"], dims="source"),
"label": xr.DataArray(["Janschwalde"], dims="source"),
"lon": xr.DataArray([14.457959], dims="source"),
"lat": xr.DataArray([51.836241], dims="source"),
"diameter": xr.DataArray([1000.0], dims="source"),
"NOx_emissions": xr.DataArray([0.42778918548939077], dims="source"), # Bottom-up estimate of emissions (in kg/s)
})
# Setup domain for plotting
NAME = "6298" # CORSO ID for Janschwalde power plant
SOURCE = sources.sel(source=[NAME])
DOMAIN = ddeq.misc.Domain.around_source(SOURCE)
CRS = ddeq.misc.get_opt_crs(DOMAIN)
# CSF settings going 30 km downstream and 80 km across
DMAX = 30e3
PLUME_WIDTH = 80e3
Read TROPOMI Level-2 data¶
The following code reads the TROPOMI Level-2 file (S5P_RPRO_L2__NO2____20210308T110658_20210308T124828_17622_03_020400_20221107T150153.nc"). You can download TROPOMI data following the example_download_tropomi.ipynb notebook in the ddeq library.
The ddeq.sats.read_S5P function will collect the different variables from the various groups to create a xarray.Dataset ready for ddeq:
filename = "/scratch/koer/TROPOMI/NO2/S5P_RPRO_L2__NO2____20210308T110658_20210308T124828_17622_03_020400_20221107T150153.nc"
data = ddeq.sats.read_S5P(filename, gas="NO2")
data
<xarray.Dataset> Size: 399MB
Dimensions: (scanline: 4173, ground_pixel: 450, corner: 4,
layer: 34, vertices: 2)
Coordinates:
* scanline (scanline) float64 33kB 0.0 1.0 ... 4.172e+03
* ground_pixel (ground_pixel) float64 4kB 0.0 1.0 ... 449.0
* corner (corner) float64 32B 0.0 1.0 2.0 3.0
* layer (layer) float64 272B 0.0 1.0 2.0 ... 32.0 33.0
* vertices (vertices) float64 16B 0.0 1.0
lat (scanline, ground_pixel) float32 8MB -73.55 ...
lon (scanline, ground_pixel) float32 8MB -134.6 ...
orbit int64 8B 17622
Data variables: (12/16)
time (scanline) datetime64[ns] 33kB 2021-03-08T11...
time_utc (scanline) <U27 451kB '2021-03-08T11:28:32.2...
qa_value (scanline, ground_pixel) float32 8MB 0.0 ......
NO2 (scanline, ground_pixel) float32 8MB nan ......
NO2_std (scanline, ground_pixel) float32 8MB nan ......
averaging_kernel (scanline, ground_pixel, layer) float32 255MB ...
... ...
solar_zenith_angle (scanline, ground_pixel) float32 8MB ...
latitude_bounds (scanline, ground_pixel, corner) float32 30MB ...
longitude_bounds (scanline, ground_pixel, corner) float32 30MB ...
surface_altitude (scanline, ground_pixel) float32 8MB ...
surface_pressure (scanline, ground_pixel) float32 8MB ...
cloud_fraction_crb (scanline, ground_pixel) float32 8MB ...
Attributes:
original file name: S5P_RPRO_L2__NO2____20210308T110658_20210308T124828_...
data source: https://s5phub.copernicus.eu/dhus/
time: 2021-03-08T11:57:44.461487658Iterate over Level-2 orbit to find Janschwalde¶
The following code will iterate over small chunks of the Level-2 orbit to find the Janschwalde power plant in the swath. The plume area will be identified using the ERA-5 wind vector.
This requires ERA-5 data on pressure levels and on single level, which can be downloade following the example_tropomi_Matimba.ipynb and example_effective_winds.ipynb notebook in the ddeq library.
for chunk in ddeq.sats.iter_S5P_swath(filename, gas="NO2", preprocessed=False):
# Plume detection using wind vector
time = pd.Timestamp(chunk.attrs['time'])
lvl_filename = time.strftime('/input/ERA5/CORSO/raw/ERA5-lvl-%Y%m%dt0000.nc')
sng_filename = time.strftime('/input/ERA5/CORSO/raw/ERA5-sng-%Y%m%dt0000.nc')
winds = ddeq.era5.read(
sng_filename,
lvl_filename,
method='gnfra',
times=time,
extent=DOMAIN.extent,
sources=SOURCE
)
chunk = ddeq.dplume.detect_from_wind(chunk, SOURCE, winds, dmax=DMAX, crs=CRS)
if "x_nodes" in chunk: # center line was computed, i.e. Janschwalde is inside chunk
break
# Prepare local coordinate system, plume area and convert from µmol/m² to kg/m²
data = ddeq.curves.compute_natural_coords(chunk)
data = ddeq.curves.compute_plume_areas(data, plume_width=PLUME_WIDTH)
ddeq.emissions.convert_units(data, "NO2", "NO2")
Quantify emissions using NO2 observations¶
The first estimate of NOx emissions is obtained by applying the CSF method with standard values for NOx correction:
# Estimate emissions using CSF method (first estimate without any corrections)
emissions = ddeq.csf.estimate_emissions(
data,
winds,
SOURCE,
'NO2',
xmin=1e3,
f_model=1.32, # NO2-to-NOx conversion factor
decay_time=4.0 * 60**2, # NOX lifetime (4 hours)
crs=CRS,
variable="NO2_mass",
background="linear"
)
# Visualize the result
fig = ddeq.vis.plot_csf_result_compact(data, winds, emissions, SOURCE, NAME, sources, crs=CRS, domain=DOMAIN)
Air mass correction¶
AMF correction uses the NO$_2$ enhancement from the first Gaussian curve to update the NO$_2$ profiles from the TROPOMI AUX dataset. In this example, ddeq looks for the file "S5P_OPER_AUX_CTMANA_20210308T000000_20210309T000000_20221231T072039.nc" in the provided path.
# Air mass factor correction using enhancements from first Gaussian plume
data = ddeq.amf.correct_amf_from_csf_method(
data,
emissions,
source_name=NAME,
path="/input/CORSO/TROPOMI/aux/"
)
# 2nd emission quantification after AMF correction
emissions = ddeq.csf.estimate_emissions(
data,
winds,
SOURCE,
'NO2',
xmin=1e3,
f_model=1.32,
decay_time=4.0 * 60**2,
crs=CRS,
variable="NO2_amf_mass",
background="linear"
)
fig = ddeq.vis.plot_csf_result_compact(data, winds, emissions, SOURCE, NAME, sources, crs=CRS, domain=DOMAIN)
NOx correction and uncertainty estimation¶
The following code applies the local NOx correction using the empirical formula from Meier et al. 2024. The code also conducts a Monte Carlo simulation to estimate the uncertainty of the estimated emissions.
# Prepare emissions dataset (as code below assumes a time series)
emissions = xr.concat([emissions.copy()], dim="time")
N = emissions.time.size
# Line density ensemble
q = np.squeeze(emissions.NO2_line_density.values)
q_std = np.squeeze(emissions.NO2_line_density_precision.values)
q_sample = q + q_std * np.random.randn(10000)
# Wind speed ensemble using Log-normal distribution with mean of u_eff and RMSE of about 1.0 m/s
u = emissions.wind_speed.values.flatten()
u_std = 1.0
mu = np.log(u**2 / np.sqrt(u**2 + u_std**2))
sigma = np.sqrt(np.log(1 + u_std**2 / u**2))
u_sample = np.exp(mu + sigma * np.random.randn(10000, N))
# Compute NO2-to-NOx conversion (between 1 km and 30 km downstream of source
along = np.linspace(1, 30, 200) * 1e3
# NO2-to-NOx conversion ensemble
seconds_since_emission = along[None,:] / u[:,None]
seconds_sample = along[None,:,None] / u_sample[:,None,:]
NO2toNOx_model = ddeq.functions.NO2toNOxConversion(name="Janschwalde")
f_local, f_local_std = NO2toNOx_model(seconds_since_emission)
f_local_sample, _ = NO2toNOx_model(seconds_sample)
f_local_sample = f_local_sample.mean(1) + f_local_std.mean(1)[None,:] * np.random.randn(10000,N)
# NOx lifetime ensemble with median of 2.5 hours
tau = 2.5 * 60**2 # in seconds
tau_std = 0.8 * 60**2 # in seconds
mu = np.log(tau**2 / np.sqrt(tau**2 + tau_std**2))
sigma = np.sqrt(np.log(1 + tau_std**2 / tau**2))
tau_local_sample = np.exp(mu + sigma * np.random.randn(10000, N))
# Decay correction term
D_local_sample = np.exp(-along[None,:,None] / (tau_local_sample[:,None,:] * u_sample[:,None,:])).mean(1)
# NOx correction term
c_local_sample = f_local_sample / D_local_sample
# Create xarray dataset
sample = xr.Dataset()
sample["time"] = xr.DataArray(emissions.time.values, dims="time")
sample["model"] = xr.DataArray(["local"], dims="model")
sample["wind_speed"] = xr.DataArray(u_sample, dims=("n", "time"))
dims = ("model", "n", "time")
sample["NOx_NO2_ratio"] = xr.DataArray([f_local_sample], dims=dims, attrs={"models": ("Janschwalde",)})
sample["NOx_lifetime"] = xr.DataArray([tau_local_sample], dims=dims, attrs={"units": "s"})
sample["NOx_decay_correction"] = xr.DataArray([D_local_sample], dims=dims)
sample["NOx_correction"] = xr.DataArray([c_local_sample], dims=dims)
# Compute NOx emissions from NOx correction (c), line density (q) and wind speed (u)
u_sample = sample["wind_speed"].values
c_sample = sample["NOx_correction"].sel(model="local").values
Q_sample = c_sample * q_sample * u_sample
# Median and +/- sigma uncertainty
Q_low, Q_mid, Q_high = scipy.stats.scoreatpercentile(Q_sample, [50-34, 50, 50+34] )
f_low, f_mid, f_high = scipy.stats.scoreatpercentile(sample["NOx_NO2_ratio"], [50-34, 50, 50+34] )
# Apply update to results dataset
emissions.f[:] = f_mid
emissions.f_precision[:] = (f_high - f_low) / 2
emissions.NOx_decay_time[:] = 2.5 * 60**2
if np.squeeze(emissions.NO2_shift) > 10e3 or np.squeeze(emissions.NO2_standard_width) > 10e3:
emissions.NOx_emissions[:] = np.nan
emissions.NOx_emissions_precision[:] = np.nan
else:
emissions.NOx_emissions[:] = Q_mid
emissions.NOx_emissions_precision[:] = (Q_high - Q_low) / 2
# Plot final result
fig = ddeq.vis.plot_csf_result_compact(data, winds, emissions.isel(time=0), SOURCE, NAME, sources, crs=CRS, domain=DOMAIN)